train_grid_layer Subroutine

public subroutine train_grid_layer(kohonen_map, input_data)

Subroutine to train the grid layer of a two_level self_organized_map

Type Bound

two_level_self_organizing_map

Arguments

Type IntentOptional Attributes Name
class(two_level_self_organizing_map) :: kohonen_map

A two_level_self_organizing_map object

type(kohonen_pattern), intent(inout), dimension(:) :: input_data

A kohonen_pattern array


Calls

proc~~train_grid_layer~~CallsGraph proc~train_grid_layer two_level_self_organizing_map%train_grid_layer float float proc~train_grid_layer->float none~distance~8 kohonen_prototype%distance proc~train_grid_layer->none~distance~8 none~get_prototype kohonen_prototype%get_prototype proc~train_grid_layer->none~get_prototype none~set_prototype kohonen_prototype%set_prototype proc~train_grid_layer->none~set_prototype proc~calculate_distance_between_prototypes~2 two_level_self_organizing_map%calculate_distance_between_prototypes proc~train_grid_layer->proc~calculate_distance_between_prototypes~2 proc~calculate_u_matrix~2 two_level_self_organizing_map%calculate_u_matrix proc~train_grid_layer->proc~calculate_u_matrix~2 proc~index2position~2 index2position proc~train_grid_layer->proc~index2position~2 proc~kohonen_pattern_accessor kohonen_pattern%kohonen_pattern_accessor proc~train_grid_layer->proc~kohonen_pattern_accessor proc~position2index~2 position2index proc~train_grid_layer->proc~position2index~2 none~distance~8->none~get_prototype calculate calculate none~distance~8->calculate proc~calculate_distance_between_prototypes~2->proc~position2index~2 proc~calculate_u_matrix~2->none~distance~8

Called by

proc~~train_grid_layer~~CalledByGraph proc~train_grid_layer two_level_self_organizing_map%train_grid_layer proc~train_2lsom two_level_self_organizing_map%train_2lsom proc~train_2lsom->proc~train_grid_layer proc~train_two_level_som train_two_level_som proc~train_two_level_som->proc~train_2lsom

Variables

Type Visibility Attributes Name Initial
integer, public :: iteration
integer, public :: iepoch
integer, public :: ipattern
integer, public :: ix
integer, public :: iy
integer, public :: iz
integer, public :: jhit
integer, public :: ihit
integer, public :: khit
integer, public :: ineigh
integer, public :: jneigh
integer, public :: kneigh
integer, public :: idbg
integer, public :: number_variables
integer, public :: max_pattern
integer, public :: ierr
integer, public :: nx
integer, public :: ny
integer, public :: nz
integer, public :: ic
integer, public :: current_pos
integer, public :: cx
integer, public :: cy
integer, public :: cz
integer, public :: i
integer, public :: j
integer, public :: k
integer, public :: number_nodes
integer, public :: debug_option
integer, public :: pos
integer, public :: pos1
integer, public :: ix1
integer, public :: iy1
integer, public :: iz1
real(kind=wp), public :: distortion
real(kind=wp), public :: dist
real(kind=wp), public :: dist_hit
real(kind=wp), public :: maximum_radius
real(kind=wp), public :: minimum_radius
real(kind=wp), public :: current_radius
real(kind=wp), public :: alpha
real(kind=wp), public :: sigma2
real(kind=wp), public :: geometric_distance2
real(kind=wp), public :: h_neighborhood
real(kind=wp), public :: real_distance
real(kind=wp), public :: lambda
real(kind=wp), public :: current_distance
real(kind=wp), public :: d1
real(kind=wp), public :: d2
real(kind=wp), public :: d3
real(kind=wp), public :: term3
real(kind=wp), public :: distance_ratio
type(kohonen_prototype), public :: current_prototype
type(kohonen_prototype), public :: current_prototype1
real(kind=wp), public, dimension(kohonen_map%number_variables1, kohonen_map%number_variables2) :: current_values
real(kind=wp), public, dimension(kohonen_map%number_variables1, kohonen_map%number_variables2) :: prototype_values
real(kind=wp), public, dimension(kohonen_map%number_variables1, kohonen_map%number_variables2) :: winner_values
real(kind=wp), public, dimension(kohonen_map%number_variables1, kohonen_map%number_variables2) :: term1
real(kind=wp), public, dimension(kohonen_map%number_variables1, kohonen_map%number_variables2) :: term2
integer, public, allocatable :: pattern_index(:,:,:,:)
logical, public :: testop
integer, public :: unit_out

Source Code

   subroutine train_grid_layer(kohonen_map,input_data)
   !========================================================================================
!! Subroutine to train the grid layer of a two_level self_organized_map
      class(two_level_self_organizing_map) :: kohonen_map
!! A `two_level_self_organizing_map` object
      type(kohonen_pattern),dimension(:),intent(inout) :: input_data
!! A `kohonen_pattern` array
      integer :: iteration,iepoch,ipattern,ix,iy,iz,jhit,ihit,khit,ineigh,jneigh,kneigh
      integer :: idbg,number_variables,max_pattern,ierr,nx,ny,nz,ic,current_pos
      integer :: cx,cy,cz,i,j,k,number_nodes,debug_option,pos,pos1,ix1,iy1,iz1
      real(kind=wp) :: distortion,dist,dist_hit,maximum_radius,minimum_radius,current_radius,alpha
      real(kind=wp) :: sigma2,geometric_distance2,h_neighborhood,real_distance,lambda,current_distance,d1,d2,d3
      real(kind=wp) :: term3,distance_ratio
      type(kohonen_prototype) :: current_prototype,current_prototype1
      real(kind=wp),dimension(kohonen_map%number_variables1,&
                  kohonen_map%number_variables2) :: current_values,prototype_values,winner_values,term1,term2
      integer, allocatable :: pattern_index(:,:,:,:)
      logical :: testop
      integer :: unit_out
         
         nx=kohonen_map%parameters(1)%number_nodes_nx;
         ny=kohonen_map%parameters(1)%number_nodes_ny;
         nz=kohonen_map%parameters(1)%number_nodes_nz;
         idbg=kohonen_map%parameters(1)%idbg;
         debug_option=kohonen_map%parameters(1)%debug_level;
         if(debug_option .gt. 0) then
         open(idbg,file=trim(kohonen_map%parameters(1)%debug_file),status='unknown');
         endif
         iteration = 0
         distortion = 0.0_wp
         number_variables=kohonen_map%number_variables; 
         maximum_radius=real(max(kohonen_map%parameters(1)%number_nodes_nx,kohonen_map%parameters(1)%number_nodes_ny,&
                                 kohonen_map%parameters(1)%number_nodes_nz));
         lambda=2.0_wp*(1.0_wp/maximum_radius);
         minimum_radius=1.0_wp;
         if(kohonen_map%parameters(1)%view_flag) then
         write(*,*) 'TWO LEVEL SOM - Grid Layer: Training starting...'
         endif
         do iepoch = 1,kohonen_map%parameters(1)%number_epochs
            if(kohonen_map%parameters(1)%view_flag) then
               write(6,*) ' Starting epoch -- distortion',iepoch,' -- ',distortion
            endif
            distortion = 0.0_wp
            do ipattern = 1, kohonen_map%parameters(1)%number_patterns
               iteration = iteration + 1
   !          det best match grid unit
   
               ihit = 0;
               jhit = 0;
               khit = 0;
               dist_hit = 1.0e7
               call input_data(ipattern)%get(current_prototype);
               call current_prototype%get_prototype(current_values);
               !call kohonen_map%find_best_match_unit(current_prototype,ihit,jhit,khit,dist_hit);
               !call current_prototype%print(unit_out);
               !write(*,*) 'check= ',ipattern,ihit,jhit,khit,dist_hit
               if(debug_option .gt. 0) then
               write(idbg,*) 'Epoch,Current Pattern',iepoch,ipattern;
               call current_prototype%print(idbg);
               endif
               !
               do iz = 1, size(kohonen_map%grid,3)  
                  do iy = 1, size(kohonen_map%grid,2)
                     do ix=1, size(kohonen_map%grid,1)
                        dist = 0.0_wp
                        ! calculate distance
                        !call kohonen_map%grid(ix,iy,iz)%print();
                        dist=kohonen_map%grid(ix,iy,iz)%distance(current_prototype,&
                           kohonen_map%distance_function);
                        !write(*,*) ix,iy,iz,dist
                        !stop
                     if(debug_option .gt. 0) then
                        call kohonen_map%grid(ix,iy,iz)%print(idbg);
                        write(idbg,*) ix,iy,iz,dist
                     endif
                     dist = dist/float(number_variables);
                     if (dist .lt. dist_hit) then
                        dist_hit = dist
                        ihit = ix;
                        jhit = iy;
                        khit = iz;
                     endif
                  enddo!iz
               enddo!iy
            enddo!ix
            !write(*,*) 'epoch,ipat,i,j,k,d= ',iepoch,ipattern,ihit,jhit,dist_hit
   !            if(kohonen_map%parameters(1)%ireal == 7) then
   !               write(unit_out,*) 'current_prototype'
   !               call current_prototype%print(unit_out)
   !               write(unit_out,*) 'bmu ',ipattern,ihit,jhit,dist_hit
   !               call kohonen_map%grid(ihit,jhit,khit)%print(unit_out)           
   !            endif
   !           kohonen_map%number_patterns(ihit,jhit)=kohonen_map%number_patterns(ihit,jhit)+1;
               distortion = distortion + dist_hit;
            ! define radius
            current_radius = max(maximum_radius*real(1001-iteration)/1000.0 + 0.9999999999,1.0_wp);
            ! define learning rate
            alpha = max(kohonen_map%parameters(1)%learning_rate*(1.0_wp-real(iteration)/1000.0),0.01_wp);
            sigma2=current_radius**2
            !max(0.2*maximum_radius*(1.0_wp-real(iteration)/1000.0),1.0_wp);
            ! update prototypes
            if(debug_option .gt. 0) then
               write(idbg,*) 'Neighborhood,alpha= ',alpha
            endif
   !
   !           call kohonen_map%update_weights(current_values,ihit,jhit,khit,maximum_radius,iteration)
   !
            do ic=1,size(kohonen_map%coordinates,1)
               current_pos=position2index(ihit,jhit,khit,nx,ny);
               current_distance=kohonen_map%cells_distances(current_pos,ic)
               if(current_distance .lt. current_radius) then
                  geometric_distance2=current_distance**2;
                  call index2position(ic,nx,ny,nz,ineigh,jneigh,kneigh);
                  !write(*,*) ic,ineigh,jneigh,kneigh,ihit,jhit,khit
                  select case(trim(kohonen_map%parameters(1)%neighborhood_type))
                     case('gaussian')
                        h_neighborhood=alpha*exp(-0.5*geometric_distance2/sigma2);
                     case('bubble')
                        h_neighborhood=alpha;
                  end select
                  if(debug_option .gt. 0) then
                        write(idbg,*) ihit,jhit,khit,ineigh,jneigh,kneigh
                  endif
                  select case(trim(kohonen_map%parameters(1)%som_type));
                     case('normal_som');
                        call kohonen_map%grid(ineigh,jneigh,kneigh)%get_prototype(prototype_values);
                        prototype_values=prototype_values+h_neighborhood*(current_values-prototype_values);
                        call kohonen_map%grid(ineigh,jneigh,kneigh)%set_prototype(prototype_values);
                     case('visom');
                        call kohonen_map%grid(ineigh,jneigh,kneigh)%get_prototype(prototype_values);
                        call kohonen_map%grid(ihit,jhit,khit)%get_prototype(winner_values);
                        real_distance=sum((winner_values-prototype_values)**2);
                        if( (ineigh .eq. ihit) .and. (jneigh .eq. jhit) .and. (kneigh .eq. khit) ) then                           
                           prototype_values=prototype_values+h_neighborhood*(current_values-prototype_values);
                        else
                           distance_ratio=dsqrt(real_distance)/(dsqrt(geometric_distance2)*lambda);
                           term1=(current_values-winner_values);
                           term2=(winner_values-prototype_values);
                           !eps=max(1.0_wp*time_factor,0.0_wp);
                           term3=1.0_wp;!((1.0_wp-eps)+eps)
                           prototype_values=prototype_values+h_neighborhood*(term1+term2*(distance_ratio-1.0_wp)*term3);
                        endif
                        call kohonen_map%grid(ineigh,jneigh,kneigh)%set_prototype(prototype_values);                         
                     end select                     
               endif
            enddo
   !
            enddo !ipattern
            !if(kohonen_map%parameters(1)%ireal == 7) stop
         enddo!iepoch
         if(kohonen_map%parameters(1)%view_flag) then
            write(*,*) 'TWO LEVEL SOM - Grid Layer: Training finished'
         endif
   ! Print prototypes
   !     print prototypes
         inquire(unit=kohonen_map%parameters(1)%iprot,opened=testop);
         if(testop) then
         do iz=1,size(kohonen_map%grid,3);
            write(kohonen_map%parameters(1)%iprot,'(A,I4)') 'Layer ',iz
            do iy=1,size(kohonen_map%grid,2);
               do ix=1,size(kohonen_map%grid,1);               
                  write(kohonen_map%parameters(1)%iprot,'(A6,1X,3I4)') 'node= ',ix,iy,iz            
                  call kohonen_map%grid(ix,iy,iz)%print(kohonen_map%parameters(1)%iprot);
               enddo
            enddo
         enddo!ix
         endif
   !     calculate and print distance matrix
         inquire(unit=kohonen_map%parameters(1)%idist,opened=testop);
         if(testop) then
         call kohonen_map%calculate_distance_between_prototypes()
   !        
   !         do ix=1,size(kohonen_map%grid,1);
   !           do iy=1,size(kohonen_map%grid,2);
   !              do iz=1,size(kohonen_map%grid,3);
   !                  current_prototype=kohonen_map%grid(ix,iy,iz);
   !                  pos=ix+(iy-1)*kohonen_map%parameters(1)%number_nodes_ny+&
   !                    (iz-1)*kohonen_map%parameters(1)%number_nodes_nx*kohonen_map%parameters(1)%number_nodes_ny;
   !                  do ix1=1,size(kohonen_map%grid,1);
   !                     do iy1=1,size(kohonen_map%grid,2);
   !                        do iz1=1,size(kohonen_map%grid,3);
   !                           pos1=ix1+(iy1-1)*kohonen_map%parameters(1)%number_nodes_ny+&
   !                              (iz1-1)*kohonen_map%parameters(1)%number_nodes_nx*&
   !                              kohonen_map%parameters(1)%number_nodes_ny;
   !                           current_prototype1=kohonen_map%grid(ix1,iy1,iz1);
   !                           kohonen_map%distance(pos,pos1)=current_prototype1%distance(current_prototype,&
   !                                                         kohonen_map%distance_function);
   !                        enddo!iz1        
   !                     enddo!iy1
   !                  enddo!ix1
   !               enddo!iz
   !            enddo!iy         
   !         enddo!ix
   ! !
   !         do ix=1,size(kohonen_map%distance,1)
   !            write(kohonen_map%parameters(1)%idist,*) (kohonen_map%distance(ix,iy),iy=1,size(kohonen_map%distance,2));
   !         enddo!ix
         endif
   !
   !     final best match
         if(kohonen_map%parameters(1)%view_flag) then
         write(6,*) 
         write(6,*) 'TWO LEVEL SOM - Grid Layer: Find Best Match Unit...'
         write(6,*)
         endif
         max_pattern=0;
         do ipattern = 1, kohonen_map%parameters(1)%number_patterns
            ihit = 0;
            jhit = 0;
            khit = 0;
            dist_hit = 100000.0;
            call input_data(ipattern)%get(current_prototype);
            !call current_prototype%get_prototype(current_values);
   !         call kohonen_map%find_best_match_unit(current_prototype,ihit,jhit,khit,dist_hit);
            !write(*,*) 'bmu=',ipattern,ihit,jhit,khit
            do iz = 1, size(kohonen_map%grid,3);  
               do iy = 1, size(kohonen_map%grid,2);
                  do ix = 1,size(kohonen_map%grid,1)
                     dist = 0.0_wp;
                     dist=kohonen_map%grid(ix,iy,iz)%distance(current_prototype,kohonen_map%distance_function);
                     dist = dist/float(number_variables);
                     if (dist .lt. dist_hit) then
                        dist_hit = dist
                        ihit = ix;
                        jhit = iy;
                        khit = iz;
                     endif
                  enddo
               enddo
            enddo
            kohonen_map%number_patterns(ihit,jhit,khit)=kohonen_map%number_patterns(ihit,jhit,khit)+1;
            if(kohonen_map%number_patterns(ihit,jhit,khit) .gt. max_pattern) then 
            max_pattern=kohonen_map%number_patterns(ihit,jhit,khit);
            endif
            kohonen_map%cells_index(ipattern,1)=ihit;
            kohonen_map%cells_index(ipattern,2)=jhit;
            kohonen_map%cells_index(ipattern,3)=khit;
            
            if(debug_option .gt. 0) then
            write(idbg,*) ipattern,ihit,jhit,khit
            endif
            inquire(unit=kohonen_map%parameters(1)%iindex,opened=testop);
            if(testop) then
            write(kohonen_map%parameters(1)%iindex,*) ipattern,ihit,jhit,khit
            endif
   !          mtchx(i) = ihit
   !          mtchy(i) = jhit
   !          ioutrep(ihit,jhit) = i
   
         enddo !ipattern
         if(kohonen_map%parameters(1)%debug_level .gt. 0) then
         close(idbg)
         endif
         allocate(pattern_index(size(kohonen_map%grid,1),size(kohonen_map%grid,2),&
                  size(kohonen_map%grid,3),&
                  max_pattern),stat=ierr);
         pattern_index=-1;
         do ipattern=1,kohonen_map%parameters(1)%number_patterns
            ix=kohonen_map%cells_index(ipattern,1);
            iy=kohonen_map%cells_index(ipattern,2);
            iz=kohonen_map%cells_index(ipattern,3);
            do i=1,max_pattern
               if(pattern_index(ix,iy,iz,i) .lt. 0) then
               pattern_index(ix,iy,iz,i)=ipattern;
               exit
               endif
            enddo
         enddo!ipattern
         inquire(unit=kohonen_map%parameters(1)%isam,opened=testop);
         if(testop) then
         do iz1=1,size(kohonen_map%grid,3);
            do iy1=1,size(kohonen_map%grid,2);
               do ix1=1,size(kohonen_map%grid,1);
                  write(kohonen_map%parameters(1)%isam,'(A,3I4)') 'Node= ',ix1,iy1,iz1
                  if(kohonen_map%number_patterns(ix1,iy1,iz1) .gt. 0) then
                     write(kohonen_map%parameters(1)%isam,'(A,100I4)') 'Sample ID= ',&
                     pattern_index(ix1,iy1,iz1,1:kohonen_map%number_patterns(ix1,iy1,iz1));
                  else
                     write(kohonen_map%parameters(1)%isam,'(A,I4)') 'Sample ID= ',0
                  endif
               enddo
            enddo
         enddo
         deallocate(pattern_index);
         endif
   !     print hit counter
      inquire(unit=kohonen_map%parameters(1)%ihit,opened=testop);
      if(testop) then
         do iz=1,size(kohonen_map%grid,3)
            write(kohonen_map%parameters(1)%ihit,'(A,I4)') 'Layer ',iz
            do ix=1,size(kohonen_map%grid,1);
               write(kohonen_map%parameters(1)%ihit,'(100I6)') (kohonen_map%number_patterns(ix,iy,iz),iy=1,size(kohonen_map%grid,2));
            enddo!ix
         enddo!iz
      endif
   
   !     calculate u_matrix
         inquire(unit=kohonen_map%parameters(1)%iumat,opened=testop);
         if(testop) then
         call kohonen_map%calculate_u_matrix();
   !       do iz=1,size(kohonen_map%grid,3);
   !          do iy=1,size(kohonen_map%grid,2);
   !             do ix=1,size(kohonen_map%grid,1);
   !                dist=0.0_wp;
   !                number_nodes=0;
   !                do k=-1,1
   !                   do j=-1,1
   !                      do i=-1,1
   !                         cx=ix+i;cy=iy+j;cz=iz+k;
   !                         if( (cx .ge. 1 .and. cx .le. size(kohonen_map%grid,1)) .and. &
   !                           (cy .ge. 1 .and. cy .le. size(kohonen_map%grid,2)) .and. &
   !                           (cz .ge. 1 .and. cz .le. size(kohonen_map%grid,3))) then
   !                               !write(*,*) ix,iy,cx,cy
   !                               dist=dist+kohonen_map%grid(ix,iy,iz)%distance(kohonen_map%grid(cx,cy,cz),&
   !                                     kohonen_map%distance_function);
   !                               number_nodes=number_nodes+1;
   !                         endif
   !                      enddo
   !                   enddo!j
   !                enddo!i
   !                kohonen_map%u_matrix(ix,iy,iz)=dist/real(number_nodes);
   !             enddo!iz
   !          enddo!iy
   !       enddo!ix
   !       do iz=1,size(kohonen_map%grid,3)
   !          write(kohonen_map%parameters(1)%iumat,'(A,I4)') 'Layer ',iz
   !          do ix=1,size(kohonen_map%grid,1)
   !             write(kohonen_map%parameters(1)%iumat,*) (kohonen_map%u_matrix(ix,iy,iz),iy=1,size(kohonen_map%grid,2)); 
   !          enddo
   !       enddo
         endif
         if(kohonen_map%parameters(1)%view_flag) then
         write(6,*) 
         write(6,*) 'TWO LEVEL SOM - Grid Layer: Find Best Match Unit finished'
         write(6,*)
         endif
         
   end subroutine train_grid_layer